1 Introduction

Since read counts are summed across cells in a pseudobulk approach, modeling continuous cell-level covariates also requires a collapsing step. Here we summarize the values of a variable from a set of cells using the mean, and store the value for each cell type. Including these variables in a regression formula uses the summarized values from the corresponding cell type.

We demonstrate this feature on a lightly modified analysis of PBMCs from 8 individuals stimulated with interferon-β (Kang, et al, 2018, Nature Biotech).

2 Standard processing

Here is the code from the main vignette:

library(dreamlet)
library(muscat)
library(ExperimentHub)
library(scater)

# Download data, specifying EH2259 for the Kang, et al study
eh <- ExperimentHub()
sce <- eh[["EH2259"]]

# only keep singlet cells with sufficient reads
sce <- sce[rowSums(counts(sce) > 0) > 0, ]
sce <- sce[, colData(sce)$multiplets == "singlet"]

# compute QC metrics
qc <- perCellQCMetrics(sce)

# remove cells with few or many detected genes
ol <- isOutlier(metric = qc$detected, nmads = 2, log = TRUE)
sce <- sce[, !ol]

# set variable indicating stimulated (stim) or control (ctrl)
sce$StimStatus <- sce$stim

In many datasets, continuous cell-level variables could be mapped reads, gene count, mitochondrial rate, etc. There are no continuous cell-level variables in this dataset, so we can simulate two from a normal distribution:

sce$value1 <- rnorm(ncol(sce))
sce$value2 <- rnorm(ncol(sce))

3 Pseudobulk

Now compute the pseudobulk using standard code:

sce$id <- paste0(sce$StimStatus, sce$ind)

# Create pseudobulk
pb <- aggregateToPseudoBulk(sce,
  assay = "counts",
  cluster_id = "cell",
  sample_id = "id",
  verbose = FALSE
)

The means per variable, cell type, and sample are stored in the pseudobulk SingleCellExperiment object:

metadata(pb)$aggr_means
## # A tibble: 128 × 5
## # Groups:   cell [8]
##    cell    id       cluster   value1   value2
##    <fct>   <fct>      <dbl>    <dbl>    <dbl>
##  1 B cells ctrl101     3.96 -0.0520   0.0107 
##  2 B cells ctrl1015    4.00 -0.00341 -0.00791
##  3 B cells ctrl1016    4     0.0394  -0.00281
##  4 B cells ctrl1039    4.04  0.0498  -0.110  
##  5 B cells ctrl107     4    -0.0517  -0.0674 
##  6 B cells ctrl1244    4     0.00521  0.0245 
##  7 B cells ctrl1256    4.01 -0.0726  -0.00362
##  8 B cells ctrl1488    4.02 -0.159   -0.0622 
##  9 B cells stim101     4.09 -0.0255  -0.0936 
## 10 B cells stim1015    4.06  0.0929   0.135  
## # ℹ 118 more rows

4 Analysis

Including these variables in a regression formula uses the summarized values from the corresponding cell type. This happens behind the scenes, so the user doesn’t need to distinguish bewteen sample-level variables stored in colData(pb) and cell-level variables stored in metadata(pb)$aggr_means.

Variance partition and hypothesis testing proceeds as ususal:

form <- ~ StimStatus + value1 + value2

# Normalize and apply voom/voomWithDreamWeights
res.proc <- processAssays(pb, form, min.count = 5)

# run variance partitioning analysis
vp.lst <- fitVarPart(res.proc, form)

# Summarize variance fractions genome-wide for each cell type
plotVarPart(vp.lst, label.angle = 60)

# Differential expression analysis within each assay
res.dl <- dreamlet(res.proc, form)

# dreamlet results include coefficients for value1 and value2
res.dl
## class: dreamletResult 
## assays(8): B cells CD14+ Monocytes ... Megakaryocytes NK cells
## Genes:
##  min: 164 
##  max: 5262 
## details(7): assay n_retain ... n_errors error_initial
## coefNames(4): (Intercept) StimStatusstim value1 value2

5 Session Info

## R version 4.4.0 alpha (2024-03-27 r86216)
## Platform: aarch64-apple-darwin20
## Running under: macOS Ventura 13.6.5
## 
## Matrix products: default
## BLAS:   /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/lib/libRblas.0.dylib 
## LAPACK: /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/lib/libRlapack.dylib;  LAPACK version 3.12.0
## 
## locale:
## [1] en_US.UTF-8/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
## 
## time zone: America/New_York
## tzcode source: internal
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] muscData_1.17.0             scater_1.32.0              
##  [3] scuttle_1.14.0              ExperimentHub_2.12.0       
##  [5] AnnotationHub_3.12.0        BiocFileCache_2.12.0       
##  [7] dbplyr_2.5.0                muscat_1.18.0              
##  [9] dreamlet_1.2.0              SingleCellExperiment_1.26.0
## [11] SummarizedExperiment_1.34.0 Biobase_2.64.0             
## [13] GenomicRanges_1.56.0        GenomeInfoDb_1.40.0        
## [15] IRanges_2.38.0              S4Vectors_0.42.0           
## [17] BiocGenerics_0.50.0         MatrixGenerics_1.16.0      
## [19] matrixStats_1.2.0           variancePartition_1.34.0   
## [21] BiocParallel_1.38.0         limma_3.60.0               
## [23] ggplot2_3.5.0               BiocStyle_2.32.0           
## 
## loaded via a namespace (and not attached):
##   [1] bitops_1.0-7              httr_1.4.7               
##   [3] RColorBrewer_1.1-3        doParallel_1.0.17        
##   [5] Rgraphviz_2.48.0          numDeriv_2016.8-1.1      
##   [7] tools_4.4.0               sctransform_0.4.1        
##   [9] backports_1.4.1           utf8_1.2.4               
##  [11] R6_2.5.1                  metafor_4.6-0            
##  [13] mgcv_1.9-1                GetoptLong_1.0.5         
##  [15] withr_3.0.0               prettyunits_1.2.0        
##  [17] gridExtra_2.3             cli_3.6.2                
##  [19] sandwich_3.1-0            labeling_0.4.3           
##  [21] sass_0.4.9                KEGGgraph_1.64.0         
##  [23] SQUAREM_2021.1            mvtnorm_1.2-4            
##  [25] blme_1.0-5                mixsqp_0.3-54            
##  [27] zenith_1.6.0              parallelly_1.37.1        
##  [29] invgamma_1.1              RSQLite_2.3.5            
##  [31] generics_0.1.3            shape_1.4.6.1            
##  [33] gtools_3.9.5              dplyr_1.1.4              
##  [35] Matrix_1.7-0              metadat_1.2-0            
##  [37] ggbeeswarm_0.7.2          fansi_1.0.6              
##  [39] abind_1.4-5               lifecycle_1.0.4          
##  [41] multcomp_1.4-25           yaml_2.3.8               
##  [43] edgeR_4.2.0               mathjaxr_1.6-0           
##  [45] gplots_3.1.3.1            SparseArray_1.4.0        
##  [47] grid_4.4.0                blob_1.2.4               
##  [49] crayon_1.5.2              lattice_0.22-6           
##  [51] beachmat_2.20.0           msigdbr_7.5.1            
##  [53] annotate_1.82.0           KEGGREST_1.44.0          
##  [55] magick_2.8.3              pillar_1.9.0             
##  [57] knitr_1.45                ComplexHeatmap_2.20.0    
##  [59] rjson_0.2.21              boot_1.3-30              
##  [61] estimability_1.5          corpcor_1.6.10           
##  [63] future.apply_1.11.2       codetools_0.2-19         
##  [65] glue_1.7.0                data.table_1.15.4        
##  [67] vctrs_0.6.5               png_0.1-8                
##  [69] Rdpack_2.6                gtable_0.3.4             
##  [71] assertthat_0.2.1          cachem_1.0.8             
##  [73] xfun_0.43                 mime_0.12                
##  [75] rbibutils_2.2.16          S4Arrays_1.4.0           
##  [77] Rfast_2.1.0               coda_0.19-4.1            
##  [79] survival_3.5-8            iterators_1.0.14         
##  [81] statmod_1.5.0             TH.data_1.1-2            
##  [83] nlme_3.1-164              pbkrtest_0.5.2           
##  [85] bit64_4.0.5               filelock_1.0.3           
##  [87] progress_1.2.3            EnvStats_2.8.1           
##  [89] bslib_0.6.2               TMB_1.9.10               
##  [91] irlba_2.3.5.1             vipor_0.4.7              
##  [93] KernSmooth_2.23-22        colorspace_2.1-0         
##  [95] rmeta_3.0                 DBI_1.2.2                
##  [97] DESeq2_1.44.0             tidyselect_1.2.1         
##  [99] emmeans_1.10.0            curl_5.2.1               
## [101] bit_4.0.5                 compiler_4.4.0           
## [103] graph_1.82.0              BiocNeighbors_1.22.0     
## [105] DelayedArray_0.30.0       bookdown_0.38            
## [107] scales_1.3.0              caTools_1.18.2           
## [109] remaCor_0.0.18            rappdirs_0.3.3           
## [111] stringr_1.5.1             digest_0.6.35            
## [113] minqa_1.2.6               rmarkdown_2.26           
## [115] aod_1.3.3                 XVector_0.44.0           
## [117] RhpcBLASctl_0.23-42       htmltools_0.5.8          
## [119] pkgconfig_2.0.3           lme4_1.1-35.2            
## [121] sparseMatrixStats_1.16.0  highr_0.10               
## [123] mashr_0.2.79              fastmap_1.1.1            
## [125] rlang_1.1.3               GlobalOptions_0.1.2      
## [127] UCSC.utils_1.0.0          DelayedMatrixStats_1.26.0
## [129] farver_2.1.1              jquerylib_0.1.4          
## [131] zoo_1.8-12                jsonlite_1.8.8           
## [133] BiocSingular_1.20.0       RCurl_1.98-1.14          
## [135] magrittr_2.0.3            GenomeInfoDbData_1.2.12  
## [137] munsell_0.5.0             Rcpp_1.0.12              
## [139] babelgene_22.9            viridis_0.6.5            
## [141] EnrichmentBrowser_2.34.0  RcppZiggurat_0.1.6       
## [143] stringi_1.8.3             zlibbioc_1.50.0          
## [145] MASS_7.3-60.2             plyr_1.8.9               
## [147] listenv_0.9.1             parallel_4.4.0           
## [149] ggrepel_0.9.5             Biostrings_2.72.0        
## [151] splines_4.4.0             hms_1.1.3                
## [153] circlize_0.4.16           locfit_1.5-9.9           
## [155] reshape2_1.4.4            ScaledMatrix_1.12.0      
## [157] BiocVersion_3.19.1        XML_3.99-0.16.1          
## [159] evaluate_0.23             RcppParallel_5.1.7       
## [161] BiocManager_1.30.22       nloptr_2.0.3             
## [163] foreach_1.5.2             tidyr_1.3.1              
## [165] purrr_1.0.2               future_1.33.2            
## [167] clue_0.3-65               scattermore_1.2          
## [169] ashr_2.2-63               rsvd_1.0.5               
## [171] broom_1.0.5               xtable_1.8-4             
## [173] fANCOVA_0.6-1             viridisLite_0.4.2        
## [175] truncnorm_1.0-9           tibble_3.2.1             
## [177] lmerTest_3.1-3            glmmTMB_1.1.9            
## [179] memoise_2.0.1             beeswarm_0.4.0           
## [181] AnnotationDbi_1.66.0      cluster_2.1.6            
## [183] globals_0.16.3            GSEABase_1.66.0